Contracting cube
Contracting cubeΒΆ
In this demo we simulate a unit cube that is fixed at \(x = 0\) and free at \(x = 1\). We use a transversally isotropic material with fiber oriented in the \(x\)-direction.
import dolfin
try:
from dolfin_adjoint import (
Constant,
DirichletBC,
Expression,
UnitCubeMesh,
interpolate,
Mesh,
)
except ImportError:
from dolfin import (
Constant,
DirichletBC,
interpolate,
Expression,
UnitCubeMesh,
Mesh,
)
import pulse
from fenics_plotly import plot
pulse.iterate.logger.setLevel(10)
# Create mesh
N = 6
mesh = UnitCubeMesh(N, N, N)
# Create subdomains
class Free(dolfin.SubDomain):
def inside(self, x, on_boundary):
return x[0] > (1.0 - dolfin.DOLFIN_EPS) and on_boundary
class Fixed(dolfin.SubDomain):
def inside(self, x, on_boundary):
return x[0] < dolfin.DOLFIN_EPS and on_boundary
# Create a facet fuction in order to mark the subdomains
ffun = dolfin.MeshFunction("size_t", mesh, 2)
ffun.set_all(0)
# Mark the first subdomain with value 1
fixed = Fixed()
fixed_marker = 1
fixed.mark(ffun, fixed_marker)
# Mark the second subdomain with value 2
free = Free()
free_marker = 2
free.mark(ffun, free_marker)
# Create a cell function (but we are not using it)
cfun = dolfin.MeshFunction("size_t", mesh, 3)
cfun.set_all(0)
# Collect the functions containing the markers
marker_functions = pulse.MarkerFunctions(ffun=ffun, cfun=cfun)
# Create mictrotructure
V_f = dolfin.VectorFunctionSpace(mesh, "CG", 1)
# Fibers
f0 = interpolate(Expression(("1.0", "0.0", "0.0"), degree=1), V_f)
# Sheets
s0 = interpolate(Expression(("0.0", "1.0", "0.0"), degree=1), V_f)
# Fiber-sheet normal
n0 = interpolate(Expression(("0.0", "0.0", "1.0"), degree=1), V_f)
# Collect the mictrotructure
microstructure = pulse.Microstructure(f0=f0, s0=s0, n0=n0)
# Create the geometry
geometry = pulse.Geometry(
mesh=mesh,
marker_functions=marker_functions,
microstructure=microstructure,
)
# Use the default material parameters
material_parameters = pulse.HolzapfelOgden.default_parameters()
# Select model for active contraction
active_model = pulse.ActiveModels.active_strain
# active_model = "active_stress"
# Set the activation
activation = Constant(0.0)
# Create material
material = pulse.HolzapfelOgden(
active_model=active_model,
parameters=material_parameters,
activation=activation,
)
# Make Dirichlet boundary conditions
def dirichlet_bc(W):
V = W if W.sub(0).num_sub_spaces() == 0 else W.sub(0)
return DirichletBC(V, Constant((0.0, 0.0, 0.0)), fixed)
# Make Neumann boundary conditions
neumann_bc = pulse.NeumannBC(traction=Constant(0.0), marker=free_marker)
# Collect Boundary Conditions
bcs = pulse.BoundaryConditions(dirichlet=(dirichlet_bc,), neumann=(neumann_bc,))
# Create problem
problem = pulse.MechanicsProblem(geometry, material, bcs)
# Solve problem
pulse.iterate.iterate(problem, activation, 0.1)
2022-03-22 21:37:41,215 [935] DEBUG pulse.iterate: Control: [0.0]
2022-03-22 21:37:41,217 [935] DEBUG pulse.iterate: Target: [0.1]
2022-03-22 21:37:41,220 [935] DEBUG pulse.iterate: Intial number of steps: 5 with step size 0.02
2022-03-22 21:37:41,221 [935] INFO pulse.iterate: Iterating to target control (f_21)...
2022-03-22 21:37:41,223 [935] INFO pulse.iterate: Current control: f_21 = 0.000
2022-03-22 21:37:41,223 [935] INFO pulse.iterate: Target: 0.100
2022-03-22 21:37:41,225 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:37:41,227 [935] DEBUG pulse.iterate: Maximum difference: 1.000e-01
2022-03-22 21:37:41,228 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:37:41,229 [935] DEBUG pulse.iterate: Current control: 0.020
2022-03-22 21:40:04,167 [935] DEBUG pulse.iterate:
SUCCESFULL STEP:
2022-03-22 21:40:04,169 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:04,170 [935] DEBUG pulse.iterate: Maximum difference: 8.000e-02
2022-03-22 21:40:04,172 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:04,173 [935] DEBUG pulse.iterate: Maximum difference: 8.000e-02
2022-03-22 21:40:04,175 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:40:04,177 [935] DEBUG pulse.iterate: Current control: 0.040
2022-03-22 21:40:08,908 [935] DEBUG pulse.iterate:
SUCCESFULL STEP:
2022-03-22 21:40:08,909 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:08,912 [935] DEBUG pulse.iterate: Maximum difference: 6.000e-02
2022-03-22 21:40:08,913 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:08,915 [935] DEBUG pulse.iterate: Maximum difference: 6.000e-02
2022-03-22 21:40:08,915 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:40:08,918 [935] DEBUG pulse.iterate: Current control: 0.060
2022-03-22 21:40:11,306 [935] DEBUG pulse.iterate: Solver did not converge...
2022-03-22 21:40:11,308 [935] DEBUG pulse.iterate:
NOT CONVERGING
2022-03-22 21:40:11,309 [935] DEBUG pulse.iterate: Reduce control step
2022-03-22 21:40:11,310 [935] DEBUG pulse.iterate: Assign old state
2022-03-22 21:40:11,317 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:11,318 [935] DEBUG pulse.iterate: Maximum difference: 6.000e-02
2022-03-22 21:40:11,319 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:40:11,320 [935] DEBUG pulse.iterate: Current control: 0.050
*** Warning: PETSc SNES solver diverged in 3 iterations with divergence reason DIVERGED_FNORM_NAN.
2022-03-22 21:40:14,218 [935] DEBUG pulse.iterate:
SUCCESFULL STEP:
2022-03-22 21:40:14,219 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:14,220 [935] DEBUG pulse.iterate: Maximum difference: 5.000e-02
2022-03-22 21:40:14,221 [935] DEBUG pulse.iterate: Adapt step size. New step size: 0.015
2022-03-22 21:40:14,222 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:14,223 [935] DEBUG pulse.iterate: Maximum difference: 5.000e-02
2022-03-22 21:40:14,225 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:40:14,226 [935] DEBUG pulse.iterate: Current control: 0.065
2022-03-22 21:40:16,610 [935] DEBUG pulse.iterate:
SUCCESFULL STEP:
2022-03-22 21:40:16,612 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:16,613 [935] DEBUG pulse.iterate: Maximum difference: 3.500e-02
2022-03-22 21:40:16,614 [935] DEBUG pulse.iterate: Adapt step size. New step size: 0.022
2022-03-22 21:40:16,616 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:16,617 [935] DEBUG pulse.iterate: Maximum difference: 3.500e-02
2022-03-22 21:40:16,618 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:40:16,619 [935] DEBUG pulse.iterate: Current control: 0.087
2022-03-22 21:40:18,402 [935] DEBUG pulse.iterate: Solver did not converge...
2022-03-22 21:40:18,404 [935] DEBUG pulse.iterate:
NOT CONVERGING
2022-03-22 21:40:18,405 [935] DEBUG pulse.iterate: Reduce control step
2022-03-22 21:40:18,407 [935] DEBUG pulse.iterate: Assign old state
2022-03-22 21:40:18,415 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:18,417 [935] DEBUG pulse.iterate: Maximum difference: 3.500e-02
2022-03-22 21:40:18,419 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:40:18,419 [935] DEBUG pulse.iterate: Current control: 0.076
*** Warning: PETSc SNES solver diverged in 3 iterations with divergence reason DIVERGED_DTOL.
2022-03-22 21:40:20,731 [935] DEBUG pulse.iterate:
SUCCESFULL STEP:
2022-03-22 21:40:20,733 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:20,735 [935] DEBUG pulse.iterate: Maximum difference: 2.375e-02
2022-03-22 21:40:20,735 [935] DEBUG pulse.iterate: Adapt step size. New step size: 0.017
2022-03-22 21:40:20,738 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:20,738 [935] DEBUG pulse.iterate: Maximum difference: 2.375e-02
2022-03-22 21:40:20,740 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:40:20,741 [935] DEBUG pulse.iterate: Current control: 0.093
2022-03-22 21:40:23,170 [935] DEBUG pulse.iterate:
SUCCESFULL STEP:
2022-03-22 21:40:23,171 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:23,173 [935] DEBUG pulse.iterate: Maximum difference: 6.875e-03
2022-03-22 21:40:23,174 [935] DEBUG pulse.iterate: Adapt step size. New step size: 0.025
2022-03-22 21:40:23,176 [935] DEBUG pulse.iterate: Check target reached: NO
2022-03-22 21:40:23,176 [935] DEBUG pulse.iterate: Maximum difference: 6.875e-03
2022-03-22 21:40:23,179 [935] DEBUG pulse.iterate: Change step size for final iteration
2022-03-22 21:40:23,180 [935] DEBUG pulse.iterate: Try new control
2022-03-22 21:40:23,182 [935] DEBUG pulse.iterate: Current control: 0.100
2022-03-22 21:40:24,926 [935] DEBUG pulse.iterate:
SUCCESFULL STEP:
2022-03-22 21:40:24,927 [935] DEBUG pulse.iterate: Check target reached: YES!
2022-03-22 21:40:24,929 [935] DEBUG pulse.iterate: Check target reached: YES!
([Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 47),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 150),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 223),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 313),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 354),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 428),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 469),
Coefficient(FunctionSpace(Mesh(VectorElement(FiniteElement('Lagrange', tetrahedron, 1), dim=3), 0), MixedElement(VectorElement(FiniteElement('Lagrange', tetrahedron, 2), dim=3), FiniteElement('Lagrange', tetrahedron, 1))), 502)],
[Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 45),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 148),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 221),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 311),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 352),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 426),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 467),
Coefficient(FunctionSpace(None, FiniteElement('Real', None, 0)), 500)])
# Get displacement and hydrostatic pressure
u, p = problem.state.split(deepcopy=True)
V = dolfin.VectorFunctionSpace(mesh, "CG", 1)
u_int = interpolate(u, V)
new_mesh = Mesh(mesh)
dolfin.ALE.move(new_mesh, u_int)
fig = plot(mesh, show=False)
fig.add_plot(plot(new_mesh, color="red", show=False))
fig.show()